MCExtractIsoSurface

(triv_lib/mrch_run.c:83)

Prototype:

  IPObjectStruct *MCExtractIsoSurface(char *FileName,
                                      int DataType,
                                      PointType CubeDim,
                                      int Width,
                                      int Height,
                                      int Depth,
                                      int SkipFactor,
                                      CagdRType IsoVal)


Description:

Extract a polygonal iso-surface out of volumetric data file.

Parameters:

FileName: Containing the volumetric data.
DataType: Type of scalar value in volume. Can be one of: 1 - Regular ascii (seperated by while spaces. 2 - Two bytes short integer. 3 - Four bytes long integer. 4 - One byte (char) integer. 5 - Four bytes float. 6 - Eight bytes double.
CubeDim: Width, height, and depth of a single cube, in object space coordinates.
Width: Of volumetric data set.
Height: Of volumetric data set.
Depth: Of volumetric data set.
SkipFactor: ypically 1. For 2, only every second sample is considered and for i, only every i'th sample is considered, in all axes.
IsoVal: At which to extract the iso-surface.


Returned Value:

IPObjectStruct *: A polygonal approximation of the iso-surface at IsoVal computed using the marching cubes algorithm.


See Also:

MCThresholdCube MCExtractIsoSurface2

Keywords:




MCExtractIsoSurface2

(triv_lib/mrch_run.c:200)

Prototype:

  IPObjectStruct *MCExtractIsoSurface2(TrivTVStruct *TV,
                                       int Axis,
                                       CagdBType TrivarNormals,
                                       PointType CubeDim,
                                       int SkipFactor,
                                       CagdRType IsoVal)


Description:

Extract a polygonal iso-surface out of a trivariate function.

Parameters:

TV: The trivariate to compute an iso surface approximation for,
Axis: of the trivariate to handle, 1 for X, 2 for Y, etc.
TrivarNormals: If TRUE normal are computed using gradient of the trivariate, if FALSE normal are estimated using differencing.
CubeDim: Width, height, and depth of a single cube, in object space coordinates.
SkipFactor: ypically 1. For 2, only every second sample is considered and for i, only every i'th sample is considered, in all axes.
IsoVal: At which to extract the iso-surface.


Returned Value:

IPObjectStruct *: A polygonal approximation of the iso-surface at IsoVal computed using the marching cubes algorithm.


See Also:

MCThresholdCube MCExtractIsoSurface

Keywords:




MCImprovePointOnIsoSrf

(triv_lib/mrchtriv.c:165)

Prototype:

  int MCImprovePointOnIsoSrf(PointType Pt,
                             PointType CubeDim,
                             CagdRType IsoVal,
                             CagdRType Tolerance,
                             CagdRType AllowedError)


Description:

Improves a given Pt to "sit" exactly on top of the iso surface, by along the gradient of the trivariate. Assume TV has been preprocessed by MCImprovePointOnIsoSrfPrelude.

Parameters:

Pt: Position to improve.
CubeDim: Size of a single cell in the trivariate volume.
IsoVal: Of iso surface extracted from TV that Pt is approximately on.
Tolerance: Requested accuracy.
AllowedError: Maximally allowed error to be considered valid. If zero it is ignored.


Returned Value:

int: TRUE if successful, FALSE otherwise.


See Also:

MCImprovePointOnIsoSrfPrelude MCImprovePointOnIsoSrfPostlude MCExtractIsoSurface2

Keywords:




MCThresholdCube

(triv_lib/mrchcube.c:72)

Prototype:

  MCPolygonStruct *MCThresholdCube(MCCubeCornerScalarStruct *CCS,
                                   RealType Threshold)


Description:

Given 8 cube corner values (scalars), compute the polygon(s) in this cube along the isosurface at Threshold. if CCS has gradient information, it is used to approximate normals at the vertices.
                    7            K           6
                     ***********************
                   * +                   * *
               L *   +                 *   *              Vertices 0 - 7
               *     +     I         * J   *              Edges    A - L
           4 *********************** 5     *
             *       +             *       *
             *       +             *       * G
             *       + H           *       *
             *       +             *       *
             *       +             * F     *
           E *       +       C     *       *
             *       ++++++++++++++*+++++++* 2
             *   D + 3             *     *
             *   +                 *   * B
             * +                   * *
             ***********************
            0           A           1


Parameters:

CCS: The cube's dimensions/information.
Threshold: Iso surface level.


Returned Value:

MCPolygonStruct *: List of polygons (not necessarily triangles), or possibly NULL.


See Also:

MCExtractIsoSurface

Keywords:

marching cubes


TrivBspNew

(triv_lib/triv_gen.c:96)

Prototype:

  TrivTVStruct *TrivBspTVNew(int ULength,
                             int VLength,
                             int WLength,
                             int UOrder,
                             int VOrder,
                             int WOrder,
                             CagdPointType PType)


Description:


Allocates the memory required for a new Bspline trivariate.

Parameters:

ULength: Number of control points in the U direction.
VLength: Number of control points in the V direction.
WLength: Number of control points in the W direction.
UOrder: Order of trivariate in the U direction.
VOrder: Order of trivariate in the V direction.
WOrder: Order of trivariate in the W direction.
PType: Type of control points (E2, P3, etc.).


Returned Value:

TrivTVStruct *: An uninitialized freeform trivariate Bspline.


Keywords:

trivariates allocation


TrivBspTVDegreeRaise

(triv_lib/trivrais.c:175)

Prototype:

  TrivTVStruct *TrivBspTVDegreeRaise(TrivTVStruct *TV, TrivTVDirType Dir)


Description:

Returns a new Bspline surface, identical to the original but with one degree higher, in the requested direction Dir.

Parameters:

TV: To raise it degree by one.
Dir: Direction of degree raising. Either U, V or W.


Returned Value:

TrivTVStruct *: A trivariate with one degree higher in direction Dir, representing the same geometry as TV.


Keywords:

degree raising


TrivBspTVDerive

(triv_lib/triv_der.c:154)

Prototype:

  TrivTVStruct *TrivBspTVDerive(TrivTVStruct *TV, TrivTVDirType Dir)


Description:

Given a Bspline trivariate, computes its partial derivative trivariate in direction Dir. Let old control polygon be P(i), i = 0 to k-1, and Q(i) be new one then:
Q(i) = (k - 1) * (P(i+1) - P(i)) / (Kv(i + k) - Kv(i + 1)), i = 0 to k-2.


Parameters:

TV: Trivariate to differentiate.
Dir: Direction of differentiation. Either U or V or W.


Returned Value:

TrivTVStruct *: Differentiated trivariate in direction Dir. A Bspline trivariate.


Keywords:

trivariates


TrivBspTVKnotInsertNDiff

(triv_lib/triv_ref.c:86)

Prototype:

  TrivTVStruct *TrivBspTVKnotInsertNDiff(TrivTVStruct *TV,
                                         TrivTVDirType Dir,
                                         int Replace,
                                         CagdRType *t,
                                         int n)


Description:


Given a Bspline trivariate, inserts n knots with different values as defined by t. If, however, Replace is TRUE, the knot are simply replacing the current knot vector in the prescribed direction.

Parameters:

TV: Trivariate to refine according to t in direction Dir.
Replace: If TRUE t is a knot vector exaclt in the length of the knot vector in direction Dir in TV and t simply replaces than knot vector. If FALSE, the knot vector in direction Dir in TV is refined by adding all the knots in t.
t: Knot vector to refine/replace the knot vector of TV in direction Dir.
n: Length of vector t.


Returned Value:

TrivTVStruct *: The refined trivariate. A Bspline trivariate.


Keywords:

trivariates


TrivBzrNew

(triv_lib/triv_gen.c:140)

Prototype:

  TrivTVStruct *TrivBzrTVNew(int ULength,
                             int VLength,
                             int WLength,
                             CagdPointType PType)


Description:

Allocates the memory required for a new Bezier trivariate.

Parameters:

ULength: Number of control points in the U direction.
VLength: Number of control points in the V direction.
WLength: Number of control points in the W direction.
PType: Type of control points (E2, P3, etc.).


Returned Value:

TrivTVStruct *: An uninitialized freeform trivariate Bezier.


Keywords:

trivariates allocation


TrivBzrTVDegreeRaise

(triv_lib/trivrais.c:69)

Prototype:

  TrivTVStruct *TrivBzrTVDegreeRaise(TrivTVStruct *TV, TrivTVDirType Dir)


Description:

Returns a new Bezier trivariate, identical to the original but with one degree higher, in the requested direction Dir. Let old control polygon be P(i), i = 0 to k-1, and Q(i) be new one then:
                     i            k-i
Q(0) = P(0), Q(i) = --- P(i-1) + (---) P(i), Q(k) = P(k-1).
                     k             k

This is applied to all rows/cols of the trivariate.

Parameters:

TV: To raise it degree by one.
Dir: Direction of degree raising. Either U, V or W.


Returned Value:

TrivTVStruct *: A surface with one degree higher in direction Dir, representing the same geometry as TV.


Keywords:

degree raising


TrivBzrTVDerive

(triv_lib/triv_der.c:66)

Prototype:

  TrivTVStruct *TrivBzrTVDerive(TrivTVStruct *TV, TrivTVDirType Dir)


Description:

Given a Bezier trivariate, computes its partial derivative trivariate in direction Dir. Let old control polygon be P(i), i = 0 to k-1, and Q(i) be new one then:
Q(i) = (k - 1) * (P(i+1) - P(i)), i = 0 to k-2.


Parameters:

TV: Trivariate to differentiate.
Dir: Direction of differentiation. Either U or V or W.


Returned Value:

TrivTVStruct *: Differentiated trivariate in direction Dir. A Bezier trivariate.


Keywords:

trivariates


TrivCnvrtBezier2BsplineTV

(triv_lib/triv_gen.c:513)

Prototype:

  TrivTVStruct *TrivCnvrtBezier2BsplineTV(TrivTVStruct *TV)


Description:


Converts a Bezier trivariate into a Bspline trivariate by adding two open end uniform knot vectors to it.

Parameters:

TV: A Bezier trivariate to convert to a Bspline TV.


Returned Value:

TrivTVStruct *: A Bspline trivariate representing the same geometry as the given Bezier TV.


Keywords:

conversion trivariate


TrivCoerceTVTo

(triv_lib/trivcoer.c:25)

Prototype:

  TrivTVStruct *TrivCoerceTVTo(TrivTVStruct *TV, CagdPointType PType)


Description:

Coerces a trivariate to point type PType.

Parameters:

TV: To coerce to a new point type PType.
PType: New point type for TV.


Returned Value:

TrivTVStruct *: A new trivariate with PType as its point type.


Keywords:

coercion


TrivCoverIsoSurfaceUsingStrokes

(triv_lib/mrchtriv.c:277)

Prototype:

  CagdCrvStruct *TrivCoverIsoSurfaceUsingStrokes(TrivTVStruct *TV,
                                                 int NumStrokes,
                                                 int StrokeType,
                                                 CagdPType MinMaxPwrLen,
                                                 CagdRType StepSize,
                                                 CagdRType IsoVal,
                                                 CagdVType ViewDir)


Description:

Computes a coverage of an iso surface at IsoVal of the trivariate T using curves along principal curvatures.

Parameters:

TV: To cover with strokes along principal curvatures of iso surface of value IsoVal.
NumStrokes: Number of strokes to distribute on the implicit surface.
StrokeType: 1 - Draw strokes along minimal principal curvature. 2 - Draw strokes along maximal principal curvature. 3 - Draw strokes along both principal curvatures. 4 - Draw strokes along constant X planes. 5 - Draw strokes along constant Y planes. 6 - Draw strokes along constant Z planes. 7 - Draw strokes along silhouette lines. 8 - Draw strokes orthogonal to silhouette lines. 9 - Draw strokes along both silhouette lines and lines orthogonal to silhouette lines. StrokesType >= 10 equals StrokeType < 10 but also emphasizes the silhouette areas setting longer edges along silhouettes.
MinMaxPwrLen: Arc length of each stroke (randomized in between). a triplet of the form (Min, Max, Power) that determines
the length of each stroke as: Avg = (Max + Min) / 2, Dev = (Max - Min) / 2 Length = Avg + Dev * Random(0, 1)^ Pwr
StepSize: Steps to take in the piecewise linear approximation.
IsoVal: Of iso surface of TV that coverage is to be computed for.
ViewDir: Direction of view, used for silhouette computation.


Returned Value:

int: A list of curves forming the coverage or NULL if error.


See Also:

MCmprovePointOnIsoSrf MCExtractIsoSurface2

Keywords:




TrivDbg

(triv_lib/triv_dbg.c:25)

Prototype:

  void TrivDbg(void *Obj)


Description:

Prints trivariates stderr. Should be linked to programs for debugging purposes, so trivariates may be inspected from a debugger.

Parameters:

Obj: A trivariate - to be printed to stderr.


Returned Value:

void


Keywords:

debugging


TrivDescribeError

(triv_lib/triv_err.c:66)

Prototype:

  char *TrivDescribeError(TrivFatalErrorType ErrorNum)


Description:

Returns a string describing a the given error. Errors can be raised by any member of this triv library as well as other users. Raised error will cause an invokation of TrivFatalError function which decides how to handle this error. TrivFatalError can for example, invoke this routine with the error type, print the appropriate message and quit the program.

Parameters:

ErrorNum: Type of the error that was raised.


Returned Value:

char *: A string describing the error type.


Keywords:

error handling


TrivEditSingleTVPt

(triv_lib/trivedit.c:37)

Prototype:

  TrivTVStruct *TrivEditSingleTVPt(TrivTVStruct *TV,
                                   CagdCtlPtStruct *CtlPt,
                                   int UIndex,
                                   int VIndex,
                                   int WIndex,
                                   CagdBType Write)


Description:

Provides the way to modify/get a single control point into/from a trivariate.

Parameters:

TV: Trivar to be modified/query.
CtlPt: New control point to be substituted into TV. Must carry the same PType as TV if to be written to TV.
UIndex, VIndex, WIndex: n trivar TV's control mesh to substitute/query CtlPt.
Write: If TRUE CtlPt is copied into TV, if FALSE the point is copied from TV to CtlPt.


Returned Value:

TrivCrvStruct *: If Write is TRUE, the new modified curve, if WRITE is FALSE, NULL.


Keywords:

trivar editing


TrivEvalCurvature

(triv_lib/trivcurv.c:151)

Prototype:

  CagdBType TrivEvalCurvature(CagdPType Pos,
                              CagdRType *PCurv1,
                              CagdRType *PCurv2,
                              CagdVType PDir1,
                              CagdVType PDir2)


Description:

Evaluates the principal curvatures and principal directions of the isosurface at location Pos in trivariate that was preprocessed by the TrivEvalCurvaturePredule function. Arbitrary number of invokations of this function are possible once TrivEvalCurvaturePredule is called. Also one should invoke TrivEvalCurvaturePostlude once done to release all auxiliary allocated data structures.

Parameters:

Pos: Location in the parametric space of the trivariate.
PCurv1, PCurv2: The two principal curvatures computed by this function.
PDir1, PDir2: The two principal directions computed by this function.


Returned Value:

CagdBType: TRUE if successful, FALSE otherwise.


See Also:

TrivEvalCurvaturePrelude TrivEvalCurvaturePostlude TrivEvalHessian TrivEvalGradient

Keywords:




TrivEvalGradient

(triv_lib/trivcurv.c:239)

Prototype:

  CagdBType TrivEvalGradient(CagdPType Pos, CagdVType Gradient)


Description:

Evaluates the Gradient of the isosurface at location Pos in trivariate that was preprocessed by the TrivEvalCurvaturePredule function. Arbitrary number of invokations of this function are possible once TrivEvalCurvaturePredule is called. Also one should invoke TrivEvalCurvaturePostlude once done to release all auxiliary allocated data structures.

Parameters:

Pos: Location in the parametric space of the trivariate.
Gradient: The Gradient computed by this function.


Returned Value:

CagdBType: TRUE if successful, FALSE otherwise.


See Also:

TrivEvalCurvaturePrelude TrivEvalCurvaturePostlude TrivEvalCurvature TrivEvalHessian

Keywords:




TrivEvalHessian

(triv_lib/trivcurv.c:280)

Prototype:

  CagdBType TrivEvalHessain(CagdPType Pos, CagdVType Hessian[3])


Description:

Evaluates the Hessian of the isosurface at location Pos in trivariate that was preprocessed by the TrivEvalCurvaturePredule function. Arbitrary number of invokations of this function are possible once TrivEvalCurvaturePredule is called. Also one should invoke TrivEvalCurvaturePostlude once done to release all auxiliary allocated data structures.

Parameters:

Pos: Location in the parametric space of the trivariate.
Hessian: The Hessian computed by this function.


Returned Value:

CagdBType: TRUE if successful, FALSE otherwise.


See Also:

TrivEvalCurvaturePrelude TrivEvalCurvaturePostlude TrivEvalCurvature TrivEvalHessian

Keywords:




TrivEvalTVCurvaturePostlude

(triv_lib/mrchtriv.c:76)

Prototype:

  void MCImprovePointOnIsoSrfPostlude(void)


Description:

Release all allocated auxiliary trivariate derivatives, for fast marching.

Parameters:

None


Returned Value:

void


See Also:

MCImprovePointOnIsoSrfPrelude MCImprovePointOnIsoSrf

Keywords:




TrivEvalTVCurvaturePostlude

(triv_lib/trivcurv.c:49)

Prototype:

  void TrivEvalTVCurvaturePostlude(void)


Description:

Release all allocated auxiliary trivariate derivatives, for curvature analysis.

Parameters:

None


Returned Value:

void


See Also:

TrivEvalTVCurvaturePrelude TrivEvalTVCurvature TrivEvalHessian TrivEvalGradient

Keywords:




TrivEvalTVCurvaturePrelude

(triv_lib/mrchtriv.c:112)

Prototype:

  CagdBType MCImprovePointOnIsoSrfPrelude(TrivTVStruct *TV)


Description:

Prepare the necessary derivative vector fields of TV for fast marching on the trivariate, later on.

Parameters:

TV: to process and prepare for further iso surface marching.


Returned Value:

CagdBType: TRUE if successful, FALSE otherwise.


See Also:

MCImprovePointOnIsoSrfPostlude MCImprovePointOnIsoSrf

Keywords:




TrivEvalTVCurvaturePrelude

(triv_lib/trivcurv.c:89)

Prototype:

  CagdBType TrivEvalTVCurvaturePrelude(TrivTVStruct *TV)


Description:

Prepare the necessary derivative vector fields of TV for curvature processing at prescribed locations, later on.

Parameters:

TV: to process and prepare for further curvature evaluations.


Returned Value:

CagdBType: TRUE if successful, FALSE otherwise.


See Also:

TrivEvalTVCurvaturePostlude TrivEvalTVCurvature TrivEvalHessian TrivEvalGradient

Keywords:




TrivFatalError

(triv_lib/triv_ftl.c:27)

Prototype:

  void TrivFatalError(TrivFatalErrorType ErrID)


Description:

Trap Triv_lib errors right here. Provides a default error handler for the triv library. Gets an error description using TrivDescribeError, prints it and exit the program using exit.

Parameters:

ErrID: Error type that was raised.


Returned Value:

void


Keywords:

error handling


TrivInterpTrivar

(triv_lib/trinterp.c:26)

Prototype:

  TrivTVStruct *TrivInterpTrivar(TrivTVStruct *TV)


Description:

Interpolates control points of given trivariate, preserving the order and continuity of the original trivariate.

Parameters:

TV: Trivariate to interpolate its control mesh.


Returned Value:

TrivTVStruct *: The interpolating trivariate.


Keywords:

interpolation


TrivLoadVolumeIntoT

(triv_lib/mrch_run.c:639)

Prototype:

  TrivTVStruct *TrivLoadVolumeIntoTV(char *FileName,
                                     int DataType,
                                     VectorType VolSize,
                                     VectorType Orders)


Description:

Loads a volumetric data set as a trivariate function of prescribed orders. Uniform open end conditions are created for it.

Parameters:

FileName: To load the trivariate data set from.
DataType: Type of scalar value in volume data file. Can be one of: 1 - Regular ascii (seperated by while spaces. 2 - Two bytes short integer. 3 - Four bytes long integer. 4 - One byte (char) integer. 5 - Four bytes float. 6 - Eight bytes double.
VolSize: Dimensions of trivariate voluem.
Orders: Orders of the three axis of the volume (in U, V, and W).


Returned Value:

TrivTVStrcut *: Loaded trivariate, or NULL if error.


Keywords:




TrivMakeTVsCompatible

(triv_lib/trivcmpt.c:50)

Prototype:

  CagdBType TrivMakeTVsCompatible(TrivTVStruct **TV1,
                                  TrivTVStruct **TV2,
                                  CagdBType SameUOrder,
                                  CagdBType SameVOrder,
                                  CagdBType SameWOrder,
                                  CagdBType SameUKV,
                                  CagdBType SameVKV,
                                  CagdBType SameWKV)


Description:

Given two trivariates, makes them compatible by: 1. Coercing their point type to be the same. 2. Making them have the same curve type. 3. Raising the degree of the lower one to be the same as the higher. 4. Refining them to a common knot vector (If Bspline and SameOrder). Note 3 is performed if SameOrder TRUE, 4 if SameKV TRUE. Both trivariates are modified IN PLACE.

Parameters:

TV1, TV2: Two surfaces to be made compatible, in place.
SameUOrder: If TRUE, this routine make sure they share the same U order.
SameVOrder: If TRUE, this routine make sure they share the same order.
SameWOrder: If TRUE, this routine make sure they share the same W order.
SameUKV: If TRUE, this routine make sure they share the same U knot vector and hence continuity. *
SameVKV: If TRUE, this routine make sure they share the same knot vector and hence continuity.
SameWKV: If TRUE, this routine make sure they share the same W knot vector and hence continuity.


Returned Value:

CagdBType: TRUE if successful, FALSE otherwise.


Keywords:

compatibility


TrivParamInDomain

(triv_lib/triv_aux.c:82)

Prototype:

  CagdBType TrivParamInDomain(TrivTVStruct *TV, CagdRType t, TrivTVDirType Dir)


Description:

Given a tri-variate and a domain - validate it.

Parameters:

TV: To make sure t is in its Dir domain.
t: Parameter value to verify.
Dir: Direction. Either U or V or W.


Returned Value:

CagdBType: TRUE if in domain, FALSE otherwise.


Keywords:

trivariates


TrivParamsInDomain

(triv_lib/triv_aux.c:118)

Prototype:

  CagdBType TrivParamsInDomain(TrivTVStruct *TV,
                               CagdRType u,
                               CagdRType v,
                               CagdRType w)


Description:

Given a tri-variate and a domain - validate it.

Parameters:

TV: To make sure (u, v, w) is in its domain.
u, v, w: To verify if it is in TV's parametric domain.


Returned Value:

CagdBType: TRUE if in domain, FALSE otherwise.


Keywords:

trivariates


TrivPlaneFrom4Points

(triv_lib/geomat4d.c:44)

Prototype:

  int TrivPlaneFrom4Points(TrivPType Pt1,
                           TrivPType Pt2,
                           TrivPType Pt3,
                           TrivPType Pt4,
                           TrivPlaneType Plane)


Description:

Computes a hyperplane in four space through the given four points. Based on a direct solution in Maple of:
 with(linalg);
 readlib(C);

 d := det( matrix( [ [x - x1, y - y1, z - z1, w - w1],
                     [x2 - x1, y2 - y1, z2 - z1, w2 - w1],
                     [x3 - x2, y3 - y2, z3 - z2, w3 - w2],
                     [x4 - x3, y4 - y3, z4 - z3, w4 - w3]] ) );
 coeff( d, x );
 coeff( d, y );
 coeff( d, z );
 coeff( d, w );


Parameters:

Pt1, Pt2, Pt3, Pt4: The four points the plane should go through.
Plane: Where the result should be placed.


Returned Value:

int: TRUE if successful, FALSE otherwise.


Keywords:

hyper plane plane


TrivSrfFromMesh

(triv_lib/triveval.c:381)

Prototype:

  CagdSrfStruct *TrivSrfFromMesh(TrivTVStruct *TV,
                                 int Index,
                                 TrivTVDirType Dir)


Description:


Extract a bivariate surface out of the given trivariate's mesh. The provided (zero based) Index specifies which Index to extract.

Parameters:

TV: Trivariate to extract a bivariate surface out of its mesh.
Index: Index of row/column/level of TV's mesh in direction Dir.
Dir: Direction of isosurface extraction. Either U or V or W.


Returned Value:

CagdSrfStruct *: A bivariate surface which was extracted from TV's Mesh. This surface is not necessarily on TV.


Keywords:

trivariates


TrivSrfFromTV

(triv_lib/triveval.c:202)

Prototype:

  CagdSrfStruct *TrivSrfFromTV(TrivTVStruct *TV,
                               CagdRType t,
                               TrivTVDirType Dir)


Description:

Extract an isoparametric surface out of the given tensor product trivariate. Operations should favor the CONST_W_DIR, in which the extraction is somewhat faster, if it is possible.

Parameters:

TV: To extract an isoparametric surface from at parameter value t in direction Dir.
t: Parameter value at which to extract the isosurface.
Dir: Direction of isosurface extraction. Either U or V or W.


Returned Value:

CagdSrfStruct *: A bivariate surface which is an isosurface of TV.


Keywords:

trivariates


TrivSrfToMesh

(triv_lib/triveval.c:510)

Prototype:

  void TrivSrfToMesh(CagdSrfStruct *Srf,
                     int Index,
                     TrivTVDirType Dir,
                     TrivTVStruct *TV)


Description:

Substitute a bivariate surface into a given trivariate's mesh. The provided (zero based) Index specifies which Index to extract.

Parameters:

Srf: Surface to substitute into the trivariate TV.
Index: Index of row/column/level of TV's mesh in direction Dir.
Dir: Direction of isosurface extraction. Either U or V or W.
TV: Trivariate to substitute a bivariate surface into its mesh.


Returned Value:

void


Keywords:

trivariates


TrivTV2CtrlMesh

(triv_lib/trivmesh.c:24)

Prototype:

  CagdPolylineStruct *TrivTV2CtrlMesh(TrivTVStruct *Trivar)


Description:

Extracts the control mesh of a surface as a list of polylines.

Parameters:

Srf: To extract a control mesh from.


Returned Value:

CagdPolylineStruct *: The control mesh of Srf.


Keywords:

control mesh


TrivTVBBox

(triv_lib/triv_aux.c:219)

Prototype:

  void TrivTVBBox(TrivTVStruct *Trivar, CagdBBoxStruct *BBox)


Description:

Computes a bounding box for a trivariate freeform function.

Parameters:

Trivar: To compute a bounding box for.
BBox: Where bounding information is to be saved.


Returned Value:

void


Keywords:

bbox bounding box


TrivTVCopy

(triv_lib/triv_gen.c:167)

Prototype:

  TrivTVStruct *TrivTVCopy(TrivTVStruct *TV)


Description:

Allocates and duplicates all slots of a trivariate structure.

Parameters:

TV: Trivariate to duplicate


Returned Value:

TrivTVStruct *: Duplicated trivariate.


Keywords:

trivariates


TrivTVCopyList

(triv_lib/triv_gen.c:231)

Prototype:

  TrivTVStruct *TrivTVCopyList(TrivTVStruct *TVList)


Description:

Duplicates a list of trivariate structures.

Parameters:

TVList: List of trivariates to duplicate.


Returned Value:

TrivTVStruct *: Duplicated list of trivariates.


Keywords:

trivariates


TrivTVDegreeRaise

(triv_lib/trivrais.c:35)

Prototype:

  TrivTVStruct *TrivTVDegreeRaise(TrivTVStruct *TV, TrivTVDirType Dir)


Description:

Returns a new trivariate representing the same curve as TV but with its degree raised by one.

Parameters:

TV: To raise its degree.
Dir: Direction of degree raising. Either U, V or W.


Returned Value:

TrivTVStruct *: A surface with same geometry as TV but with one degree higher.


Keywords:

degree raising


TrivTVDerive

(triv_lib/triv_der.c:35)

Prototype:

  TrivTVStruct *TrivTVDerive(TrivTVStruct *TV, TrivTVDirType Dir)


Description:

Given a trivariate, computes its partial derivative trivariate in direction Dir.

Parameters:

TV: Trivariate to differentiate.
Dir: Direction of differentiation. Either U or V or W.


Returned Value:

TrivTVStruct *: Differentiated trivariate in direction Dir.


Keywords:

trivariates


TrivTVDomain

(triv_lib/triv_aux.c:33)

Prototype:

  void TrivTVDomain(TrivTVStruct *TV,
                    CagdRType *UMin,
                    CagdRType *UMax,
                    CagdRType *VMin,
                    CagdRType *VMax,
                    CagdRType *WMin,
                    CagdRType *WMax)


Description:

Given a tri-variate, returns its parametric domain.

Parameters:

TV: Trivariate function to consider.
UMin, UMax: U Domain of TV will be placed herein.
VMin, VMax: V Domain of TV will be placed herein.
WMin, WMax: W Domain of TV will be placed herein.


Returned Value:

void


Keywords:

trivariates


TrivTVEval

(triv_lib/triveval.c:44)

Prototype:

  CagdRType *TrivTVEval(TrivTVStruct *TV, CagdRType u, CagdRType v, CagdRType w)


Description:

Evaluates the given tensor product trivariate at a given point, by extracting an isoparamteric surface along w and evaluating (u,v) in it.
          +-----------------------+
      W  /                       /|
     /  /                       / |
    /  /      U -->            /  |  Orientation of
      +-----------------------+   |  control tri-variate mesh.
  V | |P0                 Pi-1|   +
    v |Pi                P2i-1|  /
      |                       | /
      |Pn-i               Pn-1|/
      +-----------------------+


Parameters:

TV: To evaluate at given (u, v, w) parametric location.
u, v, w: Parametric location to evaluate TV at.


Returned Value:

CagdRType *: A vector holding all the coefficients of all components of the trivariate's point type. If for example trivariate point type is P2, the W, X, and Y will be saved in the first three locations of the returned vector. The first location (index 0) of the returned vector is reserved for the rational coefficient W and XYZ always starts at second location of the returned vector (index 1).


Keywords:

evaluation trivariates


TrivTVEval2

(triv_lib/triveval.c:163)

Prototype:

  CagdRType *TrivTVEval2(TrivTVStruct *TV, CagdRType u, CagdRType v, CagdRType w)


Description:


This function is the same as TrivTVEval2 above. Cleaner, but much less efficient. Evaluates the given tensor product trivariate at a given point, by extracting an isoparamteric surface along w and evaluating (u, v) in it.

Parameters:

TV: To evaluate at given (u, v, w) parametric location.
u, v, w: Parametric location to evaluate TV at.


Returned Value:

CagdRType *: A vector holding all the coefficients of all components of the trivariate's point type. If for example trivariate point type is P2, the W, X, and Y will be saved in the first three locations of the returned vector. The first location (index 0) of the returned vector is reserved for the rational coefficient W and XYZ always starts at second location of the returned vector (index 1).


Keywords:

evaluation trivariates


TrivTVFree

(triv_lib/triv_gen.c:260)

Prototype:

  void TrivTVFree(TrivTVStruct *TV)


Description:

Deallocates and frees all slots of a trivariate structure.

Parameters:

TV: Trivariate to free.


Returned Value:

void


Keywords:

trivariates


TrivTVFreeList

(triv_lib/triv_gen.c:296)

Prototype:

  void TrivTVFreeList(TrivTVStruct *TVList)


Description:

Deallocates and frees a list of trivariate structures.

Parameters:

TVList: Trivariate list to free.


Returned Value:

void


Keywords:

trivariates


TrivTVFromSrfs

(triv_lib/trivstrv.c:33)

Prototype:

  TrivTVStruct *TrivTVFromSrfs(CagdSrfStruct *SrfList, int OtherOrder)


Description:

Constructs a trivariate using a set of surfaces. Surfaces are made to be compatible and then each is substituted into the new trivariate's mesh as a row. If the OtherOrder is less than the number of curves, number of curves is used. A knot vector is formed with uniform open end for the other direction, so it interpolates the first and last surfaces. Note, however, that only the first and the last surfaces are interpolated if OtherOrder is greater than 2.

Parameters:

SrfList: List of surfaces to consturct a trivariate with.
OtherOrder: Other, third, order of trivariate.


Returned Value:

TrivTVStruct *: Constructed trivariate from surfaces.


Keywords:

trivar constructors


TrivTVListBBox

(triv_lib/triv_aux.c:246)

Prototype:

  void TrivTVListBBox(TrivTVStruct *TVs, CagdBBoxStruct *BBox)


Description:

Computes a bounding box for a list of trivariate freeform function.

Parameters:

Trivars: To compute a bounding box for.
BBox: Where bounding information is to be saved.


Returned Value:

void


Keywords:

bbox bounding box


TrivTVMatTransform

(triv_lib/triv_gen.c:481)

Prototype:

  void TrivTVMatTransform(TrivTVStruct *TV, CagdMType Mat)


Description:

Transforms, in place, the given TV as specified by homogeneous matrix Mat.

Parameters:

TV: Trivariate to transform.
Mat: Homogeneous transformation to apply to TV.


Returned Value:

void


Keywords:

trivariates


TrivTVNew

(triv_lib/triv_gen.c:35)

Prototype:

  TrivTVStruct *TrivTVNew(TrivGeomType GType,
                          CagdPointType PType,
                          int ULength,
                          int VLength,
                          int WLength)


Description:

Allocates the memory required for a new trivariate.

Parameters:

GType: Type of geometry the curve should be - Bspline, Bezier etc.
PType: Type of control points (E2, P3, etc.).
ULength: Number of control points in the U direction.
VLength: Number of control points in the V direction.
WLength: Number of control points in the W direction.


Returned Value:

TrivTVStruct *: An uninitialized freeform trivariate.


Keywords:

trivariates allocation


TrivTVRefineAtParams

(triv_lib/triv_ref.c:41)

Prototype:

  TrivTVStruct *TrivTVRefineAtParams(TrivTVStruct *TV,
                                     TrivTVDirType Dir,
                                     CagdBType Replace,
                                     CagdRType *t,
                                     int n)


Description:

Given a trivariate, refines it at the given n knots as defined by the vector t. If Replace is TRUE, the values replace the current knot vector. Returns pointer to refined TV (Note a Bezier trivariate will be converted into a Bspline trivariate).

Parameters:

TV: Trivariate to refine according to t in direction Dir.
Dir: Direction of refinement. Either U or V or W.
Replace: If TRUE t is a knot vector exaclt in the length of the knot vector in direction Dir in TV and t simply replaces than knot vector. If FALSE, the knot vector in direction Dir in TV is refined by adding all the knots in t.
t: Knot vector to refine/replace the knot vector of TV in direction Dir.
n: Length of vector t.


Returned Value:

TrivTVStruct *: The refined trivariate. Always a Bspline trivariate.


Keywords:

trivariates


TrivTVRegionFromTV

(triv_lib/triv_aux.c:147)

Prototype:

  TrivTVStruct *TrivTVRegionFromTV(TrivTVStruct *TV,
                                   CagdRType t1,
                                   CagdRType t2,
                                   TrivTVDirType Dir)


Description:

Given a tri-variate, returns a sub-region of it.

Parameters:

TV: To extract asub-region from.
t1, t1: Domain to extract from TV, in parametric direction Dir.
Dir: Direction to extract the sub-region. Either U or V or W.


Returned Value:

TrivTVStruct *: A sub-region of TV from t1 to t2 in direction Dir.


Keywords:

trivariates


TrivTVSubdivAtParam

(triv_lib/triv_sub.c:29)

Prototype:

  TrivTVStruct *TrivTVSubdivAtParam(TrivTVStruct *TV,
                                    CagdRType t,
                                    TrivTVDirType Dir)


Description:

Given a tri-variate, subdivides it at parameter value t in direction Dir.

Parameters:

TV: Trivariate to subdivide.
t: Parameter to subdivide at.
Dir: Direction of subdivision.


Returned Value:

TrivTVStruct *: A list of two trivariates, result of the subdivision.


Keywords:

trivariates


TrivTVTransform

(triv_lib/triv_gen.c:449)

Prototype:

  void TrivTVTransform(TrivTVStruct *TV, CagdRType *Translate, CagdRType Scale)


Description:

Linearly transforms, in place, given TV as specified by Translate and Scale.

Parameters:

TV: Trivariate to transform.
Translate: Translation factor.
Scale: Scaling factor.


Returned Value:

void


Keywords:

trivariates


TrivTriangleCopy

(triv_lib/triv_gen.c:344)

Prototype:

  TrivTriangleStruct *TrivTriangleCopy(TrivTriangleStruct *Triangle)


Description:

Allocates and duplicates all slots of a triangle structure.

Parameters:

TrivTriangleStruct *: Triangle to duplicate.


Returned Value:

TrivTriangleStruct *: Duplicated triangle.


Keywords:




TrivTriangleCopyList

(triv_lib/triv_gen.c:371)

Prototype:

  TrivTriangleStruct *TrivTriangleCopyList(TrivTriangleStruct *TriangleList)


Description:

Duplicates a list of triangle structures.

Parameters:

TriangleList: List of triangle to duplicate.


Returned Value:

TrivTriangleStruct *: Duplicated list of triangle.


Keywords:




TrivTriangleFree

(triv_lib/triv_gen.c:400)

Prototype:

  void TrivTriangleFree(TrivTriangleStruct *Triangle)


Description:

Deallocates and frees all slots of a triangle structure.

Parameters:

Triangle: Triangle to free.


Returned Value:

void


Keywords:




TrivTriangleFreeList

(triv_lib/triv_gen.c:422)

Prototype:

  void TrivTriangleFreeList(TrivTriangleStruct *TriangleList)


Description:

Deallocates and frees a list of triangle structures.

Parameters:

TriangleList: Triangle list to free.


Returned Value:

void


Keywords:




TrivTriangleNew

(triv_lib/triv_gen.c:320)

Prototype:

  TrivTriangleStruct *TrivTriangleNew(void)


Description:

Allocates the memory required for a new triangle.

Parameters:

None


Returned Value:

TrivTriangleStruct *: An uninitialized triangle.


Keywords:

allocation


TrivTwoTVsMorphing

(triv_lib/trivmrph.c:36)

Prototype:

  TrivTVStruct *TrivTwoTVsMorphing(TrivTVStruct *TV1,
                                   TrivTVStruct *TV2,
                                   CagdRType Blend)


Description:

Given two compatible trivariates (See function TrivMakeTVsCompatible), computes a convex blend between them according to Blend which must be between zero and one. Returned is the new blended trivariate.

Parameters:

TV1, TV2: The two trivariates to blend.
Blend: A parameter between zero and one


Returned Value:

TrivTVStruct *: TV2 * Blend + TV1 * (1 - Blend).


See Also:

SymbTwoCrvsMorphing SymbTwoCrvsMorphingCornerCut SymbTwoCrvsMorphingMultiRes SymbTwoSrfsMorphing TrivMakeTVsCompatible

Keywords:

morphing


TrivVectCross3Vecs

(triv_lib/geomat4d.c:180)

Prototype:

  void TrivVectCross3Vecs(TrivVType A,
                          TrivVType B,
                          TrivVType C,
                          TrivVType Res)


Description:

Computes a vector in R^4 that is perpendicular to the given three vectors.
 with(linalg);
 readlib(C);

 d := det( matrix( [ [I, J, K, L],
                   [A[0], A[1], A[2], A[3]],
                   [B[0], B[1], B[2], B[3]],
                   [C[0], C[1], C[2], C[3]] ] ) );
 coeff( d, I );
 coeff( d, J );
 coeff( d, K );
 coeff( d, L );


Parameters:

A, B, C: The three vectors to compute their cross product.
Res: Where the output goes into.


Returned Value:

void


Keywords:

cross product